Modeling

Quantum model

class qtealeaves.modeling.QuantumModel(dim, lvals, name='Rydberg', input_file='model.in', map_type='HilbertCurveMap')[source]

The class represents a physical model, e.g., how to build the Hamiltonian. Therefore, multiple instances of the class _ModelTerm can be added.

Arguments

dimint

number of spatial dimensions, i.e., the information if we have a 1d, 2d, or 3d system.

lvalsint, str, callable

Information about the number of sites; at the moment regarding one spatial dimension. For example, it is the length of the side of the cube in 3d systems.

namestr, callable

(deprecated) Name tag for the model. Only use for v1 of the input processor when writing json files. It can be parameterized via a callable returning the filename or a string-key in the parameters dictionary. Default to Rydberg

input_filestr

Name of the input file within the input folder.

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations. Default to HilbertCurveMap

add_hterm(term)[source]

Add a new term to the model.

Arguments

terminstance of _ModelTerm

Represents the additional term in the model.

apply_ham_to_state(state, ops, params)[source]

Apply the Hamiltonian to a state by contracting all operators to it, without constructing the matrix.

Arguments

stateinstance of StateVector or numpy.ndarray

The input state.

opsinstance of TNOperators

The operators consisting the Hamiltonian.

paramsdict

The parameter dictionary for the simulation.

build_ham(ops, params, qubit_limit=None)[source]

Build Hamiltonian as a matrix, which only works for very small system sizes. Moreover, the underlying terms need to support building the Hamiltonian as a matrix.

Arguments

opsinstance of TNOperators

To build the Hamiltonian, we need the actual operators which are passed within this variable.

paramsdict

The parameter dictionary for the simulation.

qubit_limitint, optional

Prevent memory issues when build_ham accidentally used for many-body systems. Default to None, i.e., 16-qubit-equivalents for dense matrices.

build_ham_sparse(ops, params, qubit_limit=None)[source]

Build Hamiltonian as a sparse matrix, which only works for very small system sizes. Moreover, the underlying terms need to support building the Hamiltonian as a matrix.

Arguments

opsinstance of TNOperators

To build the Hamiltonian, we need the actual operators which are passed within this variable.

paramsdict

The parameter dictionary for the simulation.

collect_operators()[source]

The required operators must be provided through this method. It relies on the same method of each _ModelTerm.

density_matrix(ops, params, temp, return_vec=False, k_b=1, eps_p=1e-08, max_prob=None)[source]

Diagonalize a Hamiltonian and compute its density matrix at finite temperature.

Parameters

opsinstance of TNOperators

To build the Hamiltonian, we need the actual operators which are passed within this variable.

paramsdict

The parameter dictionary for the simulation.

tempfloat

Temperature.

return_vecBoolean, optional

If True, return eigenvectors instead of density matrix. Default to False.

k_bfloat, optional

Value for the Boltzmann constant.Default to 1.

eps_pfloat, optional

Value after which the finite temperature state probabilities are cut off. Default to 1e-8.

max_probint or None, optional

Maximum number of finite temperature probabilities left after truncating. If None, there is no limit on number of probabilities. Default to None.

Return

if return_vec is False:
rho2D np.ndarray

Density matrix.

if return_vec is True:
2D np.ndarray :

Array with eigenstates of a Hamiltonian as columns. Note that this array contains only the eigenstates whose corresponding finite temperature probabilities are larger than eps_p.

prob1D np.ndarray

Finite temperature probabilities.

eval_lvals(params)[source]

Evaluate the system size via a parameter dictionary and eval_numeric_params.

Arguments

paramsdict

The parameter dictionary for the simulation.

get_coupl_str(elem, params, strength=None)[source]

Evaluate the coupling strength and return a string representation.

Arguments

eleminstance of _ModelTerm

Evaluate the coupling strength for the term in the model.

paramsdict

The parameter dictionary for the simulation.

strengthNone or numeric, optional

If None, the coupling strength will be evaluated. If given, then it will be taken as is. Default to None.

get_number_of_sites(params)[source]

Return the total number of sites in the system, i.e., the product along all dimensions.

Arguments

paramsdict

The parameter dictionary for the simulation.

get_number_of_sites_xyz(params)[source]

Returns a scalar (integer), which is the dimensions in x-, y-, and z-direction. Equal dimensions along a dimensions are assumed and only the dimension in x-direction is returned.

Arguments

paramsdict

The parameter dictionary for the simulation.

parameterization_write(folder_name_input, params, idx=0)[source]

Write the parameterization of the model into a file and return the filename.

Arguments

folder_name_inputstr

Name of the input folder, where the parameterization file will be stored.

paramsdict

The parameter dictionary for the simulation.

idxint, optional

The parameterization file allows to be indexed to support multiple files. Default to 0.

timedependent_parameter_to_dict(params)[source]

Return a dictionary with all the time dependent parameters of the simulation.

Arguments

paramsdict

The parameter dictionary for the simulation.

write_coupling_file(filename, elem, params, strength=None)[source]

Open a file and write the coupling there.

Arguments

filenamestr

Write couplings to this file. Must include the full path.

eleminstance of _ModelTerm

Evaluate the coupling strength for the term in the model.

paramsdict

The parameter dictionary for the simulation.

strengthNone or numeric, optional

If None, the coupling strength will be evaluated. If given, then it will be taken as is. Default to None.

write_input(folder_name_input, params, operator_map, mpo_mode, spo_cls)[source]

Third generation of input file writing for parameterized models. The function returns the full path to the top-level file for the model.

Arguments

folder_name_inputstr

This is the input folder where the nml file can be saved.

paramsdict

Parameters for the specific simulation.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

mpo_modeint

Signaling which MPO the user is requesting from the simulation. Choosing the MPO representation inside the simulation, where the value -1 enables an auto-selection. The other values are TPO via iTPO (0), iTPO (1), iTPO with update (2, for quenches), iTPO with compression (3), and sparse MPO (10, partially enabled), indexed sparse MPO (11, partially enabled) Default to -1.

spo_cls : constructor for the sparse MPO to be built for MPO mode 10, 11

write_sparse_mpo(folder_name_input, params, operator_map, param_map, mpo_mode, spo_cls)[source]

Write the sparse MPO to a file for a fortran simulation.

Arguments

folder_name_inputstr

This is the input folder where the file can be saved.

paramsdict

The parameter dictionary for the simulation.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapOrderedDict

For mapping parameterizations to integer IDs.

mpo_modeint

Signaling which MPO the user is requesting from the simulation. Choosing the MPO representation inside the simulation, where the value -1 enables an auto-selection. The other values are TPO via iTPO (0), iTPO (1), iTPO with update (2, for quenches), iTPO with compression (3), and sparse MPO (10, partially enabled), indexed sparse MPO (11, partially enabled) Default to -1.

spo_cls : constructor for the sparse MPO to be built for MPO mode 10, 11

Local terms

class qtealeaves.modeling.LocalTerm(operator, strength=1, prefactor=1, mask=None)[source]

Local Hamiltonian terms are versatile and probably part of any model which will be implemented. For example, the external field in the quantum Ising model can be represented as a local term.

Arguments

operatorstr

String identifier for the operator. Before launching the simulation, the python API will check that the operator is defined.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

maskcallable or None, optional

The true-false-mask allows to apply the local Hamiltonians only to specific sites, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. Default to None (all sites have a local term)

Attributes

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations.

collect_operators()[source]

The required operators must be provided through this method; thus, we return the operator in the local term.

count(params)[source]

Defines length as number of terms in fortran input file, which by now depends the presence of a mask.

Arguments

paramsdictionary

Contains the simulation parameters.

get_entries(params)[source]

Return the operator and the strength of this term.

Arguments

paramsdictionary

Contains the simulation parameters.

get_fortran_str(ll, params, operator_map, param_map, dim)[source]

Get the string representation needed to write the local terms as an plocal_type for Fortran.

Arguments

llint

Number of sites along the dimensions, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

dimint

Dimensionality of the problem, e.g., a 2d system.

get_interactions(ll, params, **kwargs)[source]

Iterator returning the local terms one-by-one, e.g., to build a Hamiltonian matrix. (In that sense, the “interaction” is obviously misleading here.)

Arguments

llint

Number of sites along the dimension, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

dimint (as keyword argument!)

Dimensionality of the problem, e.g., a 2d system.

get_sparse_matrix_operators(ll, params, operator_map, param_map, sp_ops_cls, **kwargs)[source]

Construct the sparse matrix operator for this term.

Arguments

llint

System size.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

sp_ops_clscallable (e.g., constructor)

Constructor for the sparse MPO operator to be built Has input bool (is_first_site), bool (is_last_site), bool (do_vectors).

kwargskeyword arguments

Keyword arguments are passed to get_interactions

class qtealeaves.modeling.LindbladTerm(operator, strength=1, prefactor=1, mask=None)[source]

Local Lindblad operators acting at one site are defined via this term. For the arguments see See LocalTerm.check_dim.

Details

The Lindblad equation is implemented as

\[\frac{d}{dt} \rho = -i [H, \rho] + \sum \gamma (L \rho L^{\dagger} - \frac{1}{2} \{ L^{\dagger} L, \rho \})\]
property is_oqs

Status flag if term belongs to Hamiltonian or is Lindblad.

quantum_jump_apply(state, operators, params, rand_generator, **kwargs)[source]

Apply jump with this Lindblad. Contains inplace update of state.

Arguments

state_AbstractTN

Current quantum state where jump should be applied.

operatorsTNOperators

Operator dictionary of the simulation.

paramsdict

Dictionary with parameters, e.g., to extract parameters which are not in quench or to build mask.

rand_generatorrandom number generator

Needs method random(), used to decide on jump within the sites.

kwargskeyword arguments

No keyword arguments are parsed for this term.

quantum_jump_weight(state, operators, quench, time, params, **kwargs)[source]

Evaluate the unnormalized weight for a jump with this Lindblad term.

Arguments

state_AbstractTN

Current quantum state where jump should be applied.

operatorsTNOperators

Operator dictionary of the simulation.

quenchDynamicsQuench

Current quench to evaluate time-dependent couplings.

timefloat

Time of the time evolution (accumulated dt)

paramsdict

Dictionary with parameters, e.g., to extract parameters which are not in quench or to build mask.

kwargskeyword arguments

No keyword arguments are parsed for this term.

class qtealeaves.modeling.RandomizedLocalTerm(operator, coupling_entries, strength=1, prefactor=1)[source]

Randomized local Hamiltonian terms are useful to model spinglass systems where the coupling of the local term is different for each site.

Arguments

operatorstring

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

coupling_entriesnumpy ndarray of rank-1,2,3

The coupling for the different sites. These values can only be set once and cannot be time-dependent in a time-evolution. The rank depends on the usage in 1d, 2d, or 3d systems.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

count(params)[source]

Defines length as number of terms in fortran input file, which by now depends the presence of the coupling entries.

Arguments

paramsdictionary

Contains the simulation parameters.

get_fortran_str(ll, params, operator_map, param_map, dim)[source]

Get the string representation needed to write the local terms as an plocal_type for Fortran.

Arguments

llint

Number of sites along one dimension, i.e., not the total number of sites. Assuming equal number of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

dimint

Dimensionality of the problem, e.g., a 2d system.

get_interactions(ll, params, **kwargs)[source]

See LocalTerm

Two-body interaction for 1D

class qtealeaves.modeling.TwoBodyTerm1D(operators, shift, strength=1, prefactor=1, add_complex_conjg=False, has_obc=True, mask=None)[source]

The term defines an interaction between two sites of the Hilbert space. For example, the tunneling term in the Bose-Hubbard model can be represented by this term. This class represents the 1d version.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

shiftint

Defines the distance of the interaction in compliance with the systems in higher dimensions. In the end, we iterate over all sites and apply interactions to sites (x,) and (x + shift,)

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

add_complex_conjgbool, optional

(BUG ticket #1) Aims to automatically add complex conjugated terms, e.g., for the tunneling term.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True

maskcallable or None, optional

The true-false-mask allows to apply the Hamiltonians terms only to specific sites, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. The mask is applied to the site where the left operator is acting on. Default to None (all sites have the interaction)

get_entries(params)[source]

Entries based on two coordinates; one of them is the origin (0,).

Arguments

paramsdictionary

Contains the simulation parameters.

Details

To-Do: complex-conjugate terms are wrong, i.e., they work for sigma_{i}^{+} sigma_{j}^{-} or b_{i}^{dagger} b_{j}, but not for terms as b_{j} sigma_{i}^{x}!

get_fortran_str(ll, params, operator_map, param_map)[source]

See _ModelTerm.get_fortran_str_two_body().

get_interactions(ll, params, **kwargs)[source]

Iterate over all interaction. The iterator yields a dictionary with the operator keys and the generic shift based on the origin, and a list with the coordinates in the 1d system.

get_sparse_matrix_operators(ll, params, operator_map, param_map, sp_ops_cls, **kwargs)[source]

Construct the sparse matrix operator for this term.

Arguments

llint

System size.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

sp_ops_clscallable (e.g., constructor)

Constructor for the sparse MPO operator to be built Has input bool (is_first_site), bool (is_last_site), bool (do_vectors).

kwargskeyword arguments

Keyword arguments are passed to get_interactions

iter_shifts()[source]

Return all possible shifts, which is in the case of the 1d systems exactly the defined shift.

String term for 1D

class qtealeaves.modeling.StringTerm1D(operators, strength=1, prefactor=1, has_obc=True, mask=None)[source]

The string term is applied to sites in a string with length l_x in a 1d model. If symmetries are used, none of the operators is allowed to change a symmetry sector (block diagonal).

Arguments operators : list of strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

has_obcbool, optional

Defines the boundary condition. Default to True.

maskcallable or None, optional

The true-false-mask allows to apply the Hamiltonians terms only to specific plaquettes, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. The mask is applied to the site where the lower-left corner of the plaquette operator is acting on. Default to None (all sites have the interaction)

static check_dim(dim)[source]

Only available in 1d systems.

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_fortran_str(ll, params, operator_map, param_map)[source]

Get the string representation needed to write a Fortran input file.

get_interactions(ll, params, **kwargs)[source]

Description of interactions close to the TPO formulation. It works for both Periodic and Open Boundary conditions, depending on the argument has_obc.

Arguments

ll[int]

Number of sites of the 1d lattice.

paramsdictionary

Contains the simulation parameters.

class qtealeaves.modeling.TwoBodyAllToAllTerm1D(operators, coupling_matrix, strength=1, prefactor=1)[source]

Random all-to-all two-body interaction for a one-dimensional system, e.g., as in spin glasses.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

coupling_matrixnumpy ndarray of rank-2

The coupling between the different sites in the all-to-all interactions. These values can only be set once and cannot be time-dependent in a time-evolution.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

Details

The term adds the following terms: sum_{i} sum_{j>i} strength * prefactor * coupling_mat[i, j] * A_i * B_j and only the upper triangular matrix of the coupling_mat is accessed. The hermiticity of the Hamiltonian is not ensured. Terms of the form B_j A_i with j < i are also not added.

get_fortran_str(ll, params, operator_map, param_map)[source]

Get the string representation needed to write for Fortran. This method works for any two-body interaction.

Arguments

llint

Number of sites along one dimension in the system, e.g., number of sites for one side of the square in a 2d system.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

Details

We cannot use the function _ModelTerm.get_fortran_str_two_body() as we have to modify the prefactor.

get_interactions(ll, params, **kwargs)[source]

Iterate over all interaction. The iterator yields a dictionary with the operator keys and the generic shift based on the origin, and a list with the coordinates in the 1d system.

Two-body interaction for 2D

class qtealeaves.modeling.TwoBodyTerm2D(operators, shift, strength=1, prefactor=1, isotropy_xyz=True, add_complex_conjg=False, has_obc=True, mask=None)[source]

The term defines an interaction between two sites of the Hilbert space. For example, the tunneling term in the Bose-Hubbard model can be represented by this term. This class represents the 2d version.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

shiftlist of two ints

Defines the distance of the interaction. In the end, we iterate over all sites and apply interactions to sites (x, y) and (x + shift[0], y + shift[1])

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

isotropy_xyzbool, optional

If False, the defined shift will only be applied as is. If true, we permute the defined shift to cover all spatial directions. Default to True.

add_complex_conjgbool, optional

(BUG ticket #1) Aims to automatically add complex conjugated terms, e.g., for the tunneling term.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True

maskcallable or None, optional

The true-false-mask allows to apply the Hamiltonians terms only to specific sites, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. The mask is applied to the site where the left operator is acting on. Default to None (all sites have the interaction)

Attributes

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations.

static check_dim(dim)[source]

See _ModelTerm.check_dim()

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_entries(params)[source]

Entries combine the information about possible shifts based on the origin coordination (0, 0) and the shift with the information about the operators and couplings.

Arguments

paramsdictionary

Contains the simulation parameters.

Details

To-Do: adding complex conjugates for general case beyond b bdagger case.

get_fortran_str(ll, params, operator_map, param_map)[source]

See _ModelTerm.get_fortran_str_two_body().

get_interactions(ll, params, **kwargs)[source]

These interactions are closest to the TPO description iterating over specific sites within the 1d coordinates.

Arguments

lllist of ints

Number of sites along the dimensions, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

iter_shifts()[source]

Return all possible shifts, which depends on the isotropy in the 2d case.

class qtealeaves.modeling.TwoBodyTerm2DLatticeLayout(operators, distance, layout, strength=1, prefactor=1, tolerance=1e-08)[source]

The term defines an interaction between two sites of the Hilbert space, for an arbitrary lattice layout. For example, the tunneling term in the Bose-Hubbard model can be represented by this term. This class represents the 2d version.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

distancefloat

Defines the distance of the interaction.

layoutLatticeLayout, str, callable

Instance of the LatticeLayout. The LatticeLayout class stores the positions of the (nx)x(ny) grid for an arbitrary lattice layout.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

tolerancefloat, optional

Tolerance for the distance of the interaction, i.e., the absolute distance must be smaller than the tolerance to consider the sites interacting.

Attributes

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations.

get_entries(params)[source]

The function is inherited, but we overwrite it since it has no meaning in this context.

get_fortran_str(ll, params, operator_map, param_map)[source]

See _ModelTerm.get_fortran_str_two_body().

get_interactions(ll, params, **kwargs)[source]

These interactions are closest to the TPO description iterating over specific sites within the 1d coordinates.

Arguments

llint

Number of sites along the dimensions, i.e., not the total number of sites. Assuming equal number of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

iter_shifts()[source]

The function is inherited, but we overwrite it since it has no meaning in this context.

Plaquette term for 2D

class qtealeaves.modeling.PlaquetteTerm2D(operators, strength=1, prefactor=1, has_obc=True, mask=None)[source]

The plaquette term is applied to 2x2 nearest-neighbor sites in a 2d model. If symmetries are used, none of the operators is allowed to change a symmetry sector (block diagonal).

Arguments operators : list of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True If [False, True], the topology is a strip on the x axis If [True, False], the topology is a strip on the y axis If [False,False], the topology is a thorus

maskcallable or None, optional

The true-false-mask allows to apply the Hamiltonians terms only to specific plaquettes, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. The mask is applied to the site where the lower-left corner of the plaquette operator is acting on. Default to None (all sites have the interaction)

The order of the operators is for the shifts (0,0), (0,1), (1,0), (1,1)

static check_dim(dim)[source]

Only available in 2d systems.

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_fortran_str(ll, params, operator_map, param_map)[source]

Get the string representation needed to write a Fortran input file.

get_interactions(ll, params, **kwargs)[source]

Description of interactions close to the TPO formulation. It works for both Periodic and Open Boundary conditions, depending on the argument has_obc. NOTE the order of the operators is

(x2,y2) —- (x4,y4)

(x1,y1) —- (x3,y3)

Arguments

lllist of ints

Number of sites along the dimensions, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

Block term for 2D

class qtealeaves.modeling.BlockTerm2D(operators, strength=1, prefactor=1, has_obc=True, mask=None)[source]

The block term is applied to sites in a rectangle with shape l_x*l_y in a 2d model. If symmetries are used, none of the operators is allowed to change a symmetry sector (block diagonal).

Arguments operators : 2d np.array of strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True If [False, True], the topology is a strip on the x axis If [True, False], the topology is a strip on the y axis If [False,False], the topology is a thorus

maskcallable or None, optional

The true-false-mask allows to apply the Hamiltonians terms only to specific plaquettes, i.e., with true values. The function takes the dictionary with the simulation parameters as an argument. The mask is applied to the site where the lower-left corner of the plaquette operator is acting on. Default to None (all sites have the interaction)

The order of the operators is for the shifts (0,0), (0,1), (1,0), (1,1)

static check_dim(dim)[source]

Only available in 2d systems.

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_fortran_str(ll, params, operator_map, param_map)[source]

Get the string representation needed to write a Fortran input file.

Arguments

llint

Number of sites along one dimension in the system, e.g., number of sites for one side of the square in a 2d system.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

get_interactions(ll, params, **kwargs)[source]

Description of interactions close to the TPO formulation. It works for both Periodic and Open Boundary conditions, depending on the argument has_obc. NOTE the order of the operators is

(x2,y2) —- (x4,y4)

(x1,y1) —- (x3,y3)

Arguments

lllist of ints

Number of sites along the dimensions, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

Two-body interaction for 3D

class qtealeaves.modeling.TwoBodyTerm3D(operators, shift, strength=1, prefactor=1, isotropy_xyz=True, add_complex_conjg=False, has_obc=True)[source]

The term defines an interaction between two sites of the Hilbert space. For example, the tunneling term in the Bose-Hubbard model can be represented by this term. This class represents the 2d version.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

shiftlist of three ints

Defines the distance of the interaction. In the end, we iterate over all sites and apply interactions to sites (x, y, z) and (x + shift[0], y + shift[1], z + shift[2])

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

isotropy_xyzbool, optional

If False, the defined shift will only be applied as is. If true, we permute the defined shift to cover all spatial directions. Default to True.

add_complex_conjgbool, optional

(BUG ticket #1) Aims to automatically add complex conjugated terms, e.g., for the tunneling term.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True

Attributes

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations.

static check_dim(dim)[source]

See _ModelTerm.check_dim()

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_entries(params)[source]

Entries are defined based on two coordinates, the origin chosen as (0, 0, 0) and the shift. There are no site-dependent terms.

Arguments

paramsdictionary

Contains the simulation parameters.

Details

To-Do: complex-conjugate terms are wrong, i.e., they work for sigma_{i}^{+} sigma_{j}^{-} or b_{i}^{dagger} b_{j}, but not for terms as b_{j} sigma_{i}^{x}!

get_fortran_str(ll, params, operator_map, param_map)[source]

See _ModelTerm.get_fortran_str_two_body().

get_interactions(ll, params, **kwargs)[source]

These interactions are closest to the TPO description iterating over specific sites within the 1d coordinates.

Arguments

llint

Number of sites along each dimension, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

iter_shifts()[source]

Return all possible shifts, which depends on the isotropy in the 3d case.

Two-body interaction for 3D

class qtealeaves.modeling.TwoBodyTerm3D(operators, shift, strength=1, prefactor=1, isotropy_xyz=True, add_complex_conjg=False, has_obc=True)[source]

The term defines an interaction between two sites of the Hilbert space. For example, the tunneling term in the Bose-Hubbard model can be represented by this term. This class represents the 2d version.

Arguments

operatorslist of two strings

String identifier for the operators. Before launching the simulation, the python API will check that the operator is defined.

shiftlist of three ints

Defines the distance of the interaction. In the end, we iterate over all sites and apply interactions to sites (x, y, z) and (x + shift[0], y + shift[1], z + shift[2])

strengthstr, callable, numeric (optional)

Defines the coupling strength of the local terms. It can be parameterized via a value in the dictionary for the simulation or a function. Default to 1.

prefactornumeric, scalar (optional)

Scalar prefactor, e.g., in order to take into account signs etc. Default to 1.

isotropy_xyzbool, optional

If False, the defined shift will only be applied as is. If true, we permute the defined shift to cover all spatial directions. Default to True.

add_complex_conjgbool, optional

(BUG ticket #1) Aims to automatically add complex conjugated terms, e.g., for the tunneling term.

has_obcbool or list of bools, optional

Defines the boundary condition along each spatial dimension. If scalar is given, the boundary condition along each spatial dimension is assumed to be equal. Default to True

Attributes

map_typestr, optional

Selecting the mapping from a n-dimensional system to the 1d system required for the TTN simulations.

static check_dim(dim)[source]

See _ModelTerm.check_dim()

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.

get_entries(params)[source]

Entries are defined based on two coordinates, the origin chosen as (0, 0, 0) and the shift. There are no site-dependent terms.

Arguments

paramsdictionary

Contains the simulation parameters.

Details

To-Do: complex-conjugate terms are wrong, i.e., they work for sigma_{i}^{+} sigma_{j}^{-} or b_{i}^{dagger} b_{j}, but not for terms as b_{j} sigma_{i}^{x}!

get_fortran_str(ll, params, operator_map, param_map)[source]

See _ModelTerm.get_fortran_str_two_body().

get_interactions(ll, params, **kwargs)[source]

These interactions are closest to the TPO description iterating over specific sites within the 1d coordinates.

Arguments

llint

Number of sites along each dimension, i.e., not the total number of sites. Assuming list of sites along all dimension.

paramsdictionary

Contains the simulation parameters.

iter_shifts()[source]

Return all possible shifts, which depends on the isotropy in the 3d case.

Base term

class qtealeaves.modeling._ModelTerm[source]

Abstract base class for any term in a model.

static check_dim(dim)[source]

By default, do not restrict dimensionality. Overwriting this methods, allows to restrict terms to 1d, 2d, or 3d systems.

Arguments

dimint

Dimensionality of the system, e.g., 3 for a 3-d systems.

eval_strength(params)[source]

Evaluate the strength of the parameter.

Arguments

paramsdictionary

Contains the simulation parameters.

get_fortran_str_twobody(ll, params, operator_map, param_map)[source]

Get the string representation needed to write for Fortran. This method works for any two-body interaction.

Arguments

llint

Number of sites along the dimensions in the system, e.g., number of sites for both sides of the rectangle in a 2d system.

paramsdictionary

Contains the simulation parameters.

operator_mapOrderedDict

For operator string as dictionary keys, it returns the corresponding integer IDs.

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

get_interactions(ll, params, **kwargs)[source]

Iterator returning the terms one-by-one, e.g., to build a Hamiltonian matrix. Must be overwritten by inheriting class.

get_param_repr(param_map)[source]

Get the integer identifier as a string to be written into the input files for fortran.

Arguments

param_mapdict

This dictionary contains the mapping from the parameters to integer identifiers.

get_strengths()[source]

Returns an iterator over the strenghts of the term. It is just the strength in most cases, but it can vary (for example for the KrausTerm)

property is_oqs

Status flag if term belongs to Hamiltonian or is Lindblad.

static iterate_sites(ll)[source]

Iterate sites; trivial in 1d model.

Arguments

lllist of one int

Number of sites in the chain is stored in the first entry of the list.

quantum_jump_apply(state, operators, params, rand_generator, **kwargs)[source]

Apply jump with this Lindblad.

quantum_jump_weight(state, operators, quench, time, params, **kwargs)[source]

Evaluate the unnormalized weight for a jump with this Lindblad term.

class qtealeaves.modeling._ModelTerm1D[source]

Abstract base class for any term in a 1D model.

static check_dim(dim)[source]

See _ModelTerm.check_dim()

collect_operators()[source]

All the required operators are returned to ensure that they are written by fortran.